The smuthi.utility package

utility

utility.automatic_parameter_selection

Functions that assist the user in the choice of suitable numerical simulation parameters.

smuthi.utility.automatic_parameter_selection.converge_angular_resolution(simulation, detector='extinction cross section', tolerance=0.001, max_iter=30, ax=None)

Find a suitable discretization step size for the default angular arrays used for plane wave expansions.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • max_iter (int) – Break convergence loops after that number of iterations, even if no convergence has been achieved.
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

Detector value for converged settings.

smuthi.utility.automatic_parameter_selection.converge_l_max(simulation, detector='extinction cross section', tolerance=0.001, tolerance_steps=2, max_iter=100, start_from_1=True, ax=None)

Find suitable multipole cutoff degree l_max for a given particle and simulation. The routine starts with the current l_max of the particle. The value of l_max is successively incremented in a loop until the resulting relative change in the detector value is smaller than the specified tolerance. The method updates the input particle object with the l_max value for which convergence has been achieved.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object containing the particle
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • tolerance_steps (int) – Number of consecutive steps at which the tolerance must be met during multipole truncation convergence. Default: 2
  • max_iter (int) – Break convergence loop after that number of iterations, even if no convergence has been achieved.
  • start_from_1 (logical) – If true (default), start from l_max=1. Otherwise, start from the current particle l_max.
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

A 3-tuple containing
  • detector value of converged or break-off parameter settings.
  • series of lmax values
  • the detector values for the given lmax values

smuthi.utility.automatic_parameter_selection.converge_m_max(simulation, detector='extinction cross section', tolerance=0.001, target_value=None, ax=None)

Find suitable multipole cutoff order m_max for a given particle and simulation. The routine starts with the current l_max of the particle, i.e. with m_max=l_max. The value of m_max is successively decremented in a loop until the resulting relative change in the detector value is larger than the specified tolerance. The method updates the input particle object with the so determined m_max.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object containing the particle
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • max_iter (int) – Break convergence loop after that number of iterations, even if no convergence has been achieved.
  • target_value (float) – If available (typically from preceding neff selection procedure), use as target detector value
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

Detector value of converged or break-off parameter settings.

smuthi.utility.automatic_parameter_selection.converge_multipole_cutoff(simulation, detector='extinction cross section', tolerance=0.001, tolerance_steps=2, max_iter=100, current_value=None, l_max_list=None, detector_value_list=None, converge_m=True, ax=None)

Find suitable multipole cutoff degree l_max and order m_max for all particles in a given simulation object. The method updates the input simulation object with the so determined multipole truncation values.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value
  • tolerance (float) – Relative tolerance for the detector value change
  • tolerance_steps (int) – Number of consecutive steps at which the tolerance must be met during multipole truncation convergence. Default: 2
  • max_iter (int) – Break convergence loops after that number of iterations, even if no convergence has been achieved
  • current_value (float) – If specified, skip l_max run and use this value for the resulting detector value. Otherwise, start with l_max run.
  • l_max_list (list) – If current_value was specified, the l_max run is skipped. Then, this list is returned as the second item in the returned tuple.
  • detector_value_list (list) – If current_value was specified, the l_max run is skipped. Then, this list is returned as the third item in the returned tuple.
  • converge_m (logical) – If false, only converge l_max, but keep m_max=l_max. Default is true
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

A 3-tuple containing
  • detector value of converged or break-off parameter settings.
  • series of lmax values
  • the detector values for the given lmax values

smuthi.utility.automatic_parameter_selection.converge_neff_max(simulation, detector='extinction cross section', tolerance=0.001, tolerance_factor=0.1, tolerance_steps=2, max_iter=30, neff_imag=0.01, neff_resolution=0.002, neff_max_increment=0.5, neff_max_offset=0, converge_lm=True, ax=None)

Find a suitable truncation value for the multiple scattering Sommerfeld integral contour and update the simulation object accordingly.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • tolerance_factor (float) – During neff selection, a smaller tolerance should be allowed to avoid fluctuations of the order of ~tolerance which would compromise convergence. Default: 0.1
  • tolerance_steps (int) – Number of consecutive steps at which the tolerance must be met during multipole truncation convergence. Default: 2
  • max_iter (int) – Break convergence loops after that number of iterations, even if no convergence has been achieved.
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega).
  • neff_resolution (float) – Discretization of the contour (in terms of eff. refractive index).
  • neff_max_increment (float) – Increment the neff_max parameter with that step size
  • neff_max_offset (float) – Start neff_max selection from the last estimated singularity location plus this value (in terms of effective refractive index)
  • converge_lm (logical) – If set to true, update multipole truncation during each step (this takes longer time, but is necessary for critical use cases like flat particles on a substrate)
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

Detector value for converged settings.

smuthi.utility.automatic_parameter_selection.converge_neff_resolution(simulation, detector='extinction cross section', tolerance=0.001, max_iter=30, neff_imag=0.01, neff_max=None, neff_resolution=0.01, ax=None)

Find a suitable discretization step size for the multiple scattering Sommerfeld integral contour and update the simulation object accordingly.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • max_iter (int) – Break convergence loops after that number of iterations, even if no convergence has been achieved.
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega).
  • neff_max (float) – Truncation value of contour (in terms of effective refractive index).
  • neff_resolution (float) – Discretization of the contour (in terms of eff. refractive index) - start value for iteration
  • ax (np.array of AxesSubplot) – Array of AxesSubplots where to live-plot convergence output
Returns:

Detector value for converged settings.

smuthi.utility.automatic_parameter_selection.evaluate(simulation, detector)

Run a simulation and evaluate the detector. :param simulation: simulation object :type simulation: smuthi.simulation.Simulation :param detector: Specify a method that accepts a simulation as input and returns

a float. Otherwise, type “extinction cross section” to use the extinction cross section as a detector.
Returns:The detector value (float)
smuthi.utility.automatic_parameter_selection.select_numerical_parameters(simulation, detector='extinction cross section', tolerance=0.001, tolerance_factor=0.1, tolerance_steps=2, max_iter=30, neff_imag=0.01, neff_resolution=0.01, select_neff_max=True, neff_max_increment=0.5, neff_max_offset=0, neff_max=None, select_neff_resolution=True, select_angular_resolution=False, select_multipole_cutoff=True, relative_convergence=True, show_plot=True)

Trigger automatic selection routines for various numerical parameters.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object from which parameters are read and into which results are stored.
  • detector (function or string) – Function that accepts a simulation object and returns a detector value the change of which is used to define convergence. Alternatively, use “extinction cross section” (default) to have the extinction cross section as the detector value.
  • tolerance (float) – Relative tolerance for the detector value change.
  • tolerance_factor (float) – During neff selection, a smaller tolerance should be allowed to avoid fluctuations of the order of ~tolerance which would compromise convergence. Default: 0.1
  • tolerance_steps (int) – Number of consecutive steps at which the tolerance must be met during multipole truncation convergence. Default: 2
  • max_iter (int) – Break convergence loops after that number of iterations, even if no convergence has been achieved.
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega)
  • neff_resolution (float) – Discretization of the contour (in terms of eff. refractive index) - if select_neff_resolution is true, this value will be eventually overwritten. However, it is required in any case. Default: 1e-2
  • select_neff_max (logical) – If set to true (default), the Sommerfeld integral truncation parameter neff_max is determined automatically with the help of a Cauchy convergence criterion.
  • neff_max_increment (float) – Only needed if select_neff_max is true. Step size with which neff_max is incremented.
  • neff_max_offset (float) – Only needed if select_neff_max is true. Start n_eff selection from the last estimated singularity location plus this value (in terms of effective refractive index)
  • neff_max (float) – Only needed if select_neff_max is false. Truncation value of contour (in terms of effective refractive index).
  • select_neff_resolution (logical) – If set to true (default), the Sommerfeld integral discretization parameter neff_resolution is determined automatically with the help of a Cauchy convergence criterion.
  • select_angular_resolution (logical) – If set to true, the angular resolution step for the default polar and azimuthal angles is determined automatically according to a Cauchy convergenge criterion.
  • select_multipole_cutoff (logical) – If set to true (default), the multipole expansion cutoff parameters l_max and m_max are determined automatically with the help of a Cauchy convergence criterion.
  • relative_convergence (logical) – If set to true (default), the neff_max convergence and the l_max and m_max convergence routine are performed in the spirit of relative convergence, i.e., the multipole expansion convergence is checked again for each value of the Sommerfeld integral truncation. This takes more time, but is required at least in the case of flat particles near interfaces.
smuthi.utility.automatic_parameter_selection.update_contour(simulation, neff_imag=0.005, neff_max=None, neff_max_offset=0.5, neff_resolution=0.002)

Update the default k_parallel arrays in smuthi.fields with a newly constructed Sommerfeld integral contour, and set the simulation object to use the default contour for particle coupling.

Parameters:
  • simulation (smuthi.simulation.Simulation) – Simulation object
  • neff_imag (float) – Extent of the contour into the negative imaginary direction (in terms of effective refractive index, n_eff=kappa/omega).
  • neff_max (float) – Truncation value of contour (in terms of effective refractive index).
  • neff_max_offset (float) – If no value for neff_max is specified, use the last estimated singularity location plus this value (in terms of effective refractive index).
  • neff_resolution (float) – Discretization of the contour (in terms of effective refractive index).
smuthi.utility.automatic_parameter_selection.update_lmax_mmax(simulation, l_max)

Assign the same l_max and m_max = l_max to all particles in simulation

utility.cuda

smuthi.utility.cuda.enable_gpu(enable=True)

Sets the use_gpu flag to enable/disable the use of CUDA kernels.

Parameters:enable (bool) – Set use_gpu flag to this value (default=True).

utility.logging

class smuthi.utility.logging.Logger(log_filename=None, log_to_file=True, log_to_terminal=True)

Allows to prompt messages both to terminal and to log file simultaneously. It also allows to print with indentation or to temporally mute the Logger.

fileno()
flush()
write(message)
class smuthi.utility.logging.LoggerIndented(indendatation=' ')
class smuthi.utility.logging.LoggerMuted
mute_logger = <smuthi.utility.logging.Logger object>
class smuthi.utility.logging.bcolors
BOLD = '\x1b[1m'
ENDC = '\x1b[0m'
FAIL = '\x1b[91m'
HEADER = '\x1b[95m'
OKBLUE = '\x1b[94m'
OKGREEN = '\x1b[92m'
UNDERLINE = '\x1b[4m'
WARNING = '\x1b[93m'
smuthi.utility.logging.write_blue(message)
smuthi.utility.logging.write_green(message)
smuthi.utility.logging.write_header(message)
smuthi.utility.logging.write_red(message)

utility.math

This module contains several mathematical functions.

smuthi.utility.math.dx_xh(n, x)

Derivative of \(x h_n(x)\), where \(h_n(x)\) is the spherical Hankel function.

Parameters:
  • n (int) – (n>0): Order of spherical Bessel function
  • x (array, complex or float) – Argument for spherical Hankel function
Returns:

Derivative \(\partial_x(x h_n(x))\) as array.

smuthi.utility.math.dx_xj(n, x)

Derivative of \(x j_n(x)\), where \(j_n(x)\) is the spherical Bessel function.

Parameters:
  • n (int) – (n>0): Order of spherical Bessel function
  • x (array, complex or float) – Argument for spherical Bessel function
Returns:

Derivative \(\partial_x(x j_n(x))\) as array.

smuthi.utility.math.inverse_vector_rotation(r, alpha=None, beta=None, gamma=None, euler_angles=None)
smuthi.utility.math.legendre_normalized(ct, st, lmax)

Return the normalized associated Legendre function \(P_l^m(\cos\theta)\) and the angular functions \(\pi_l^m(\cos \theta)\) and \(\tau_l^m(\cos \theta)\), as defined in A. Doicu, T. Wriedt, and Y. A. Eremin: “Light Scattering by Systems of Particles”, Springer-Verlag, 2006. Two arguments (ct and st) are passed such that the function is valid for general complex arguments, while the branch cuts are defined by the user already in the definition of st.

Parameters:
  • ct (ndarray) – cosine of theta (or kz/k)
  • st (ndarray) – sine of theta (or kp/k), need to have same dimension as ct, and st**2+ct**2=1 is assumed
  • lmax (int) – maximal multipole order
Returns:

  • ndarray plm[l, m, *ct.shape] contains \(P_l^m(\cos \theta)\). The entries of the list have same dimension as ct (and st)
  • ndarray pilm[l, m, *ct.shape] contains \(\pi_l^m(\cos \theta)\).
  • ndarray taulm[l, m, *ct.shape] contains \(\tau_l^m(\cos \theta)\).

smuthi.utility.math.legendre_normalized_numbed
smuthi.utility.math.nb_wig3jj(jj_1, jj_2, jj_3, mm_1, mm_2, mm_3)
smuthi.utility.math.rotation_matrix(alpha=None, beta=None, gamma=None, euler_angles=None)
smuthi.utility.math.spherical_hankel(n, x)
smuthi.utility.math.vector_rotation(r, alpha=None, beta=None, gamma=None, euler_angles=None)
smuthi.utility.math.wigner_D(l, m, m_prime, alpha, beta, gamma, wdsympy=False)

Computation of Wigner-D-functions for the rotation of a T-matrix

Parameters:
  • l (int) – Degree \(l\) (1, …, lmax)
  • m (int) – Order \(m\) (-min(l,mmax),…,min(l,mmax))
  • m_prime (int) – Order \(m_prime\) (-min(l,mmax),…,min(l,mmax))
  • alpha (float) – First Euler angle in rad
  • beta (float) – Second Euler angle in rad
  • gamma (float) – Third Euler angle in rad
  • wdsympy (bool) – If True, Wigner-d-functions come from the sympy toolbox
Returns:

single complex value of Wigner-D-function

smuthi.utility.math.wigner_d(l, m, m_prime, beta, wdsympy=False)

Computation of Wigner-d-functions for the rotation of a T-matrix

Parameters:
  • l (int) – Degree \(l\) (1, …, lmax)
  • m (int) – Order \(m\) (-min(l,mmax),…,min(l,mmax))
  • m_prime (int) – Order \(m_prime\) (-min(l,mmax),…,min(l,mmax))
  • beta (float) – Second Euler angle in rad
  • wdsympy (bool) – If True, Wigner-d-functions come from the sympy toolbox
Returns:

real value of Wigner-d-function

utility.memoizing

Provide functionality to store intermediate results in lookup tables (memoize)

class smuthi.utility.memoizing.Memoize(fn)

To be used as a decorator for functions that are memoized.

utility.optical_constants

Provide functionality to read optical constants in format provided by refractiveindex.info website

smuthi.utility.optical_constants.read_refractive_index_from_yaml(filename, vacuum_wavelength, units='mkm', kind=1)

Read optical constants in format provided by refractiveindex.info website.

Parameters:
  • filename (str) – path and file name for yaml data downloaded from refractiveindex.info
  • vacuum_wavelength (float or np.array) – wavelengths where refractive index data is needed
  • units (str) – units for wavelength. currently, microns (‘mkm’ or ‘um’) and nanometers (‘nm’) can be selected
  • kind (int) – order of interpolation
Returns:

A pair (or np.array of pairs) of wavelength and corresponding refractive index (complex)